hysop.numerics.fft.fft module

Base interface for fast Fourier Transforms.

FFTI FFTPlanI FFTQueueI

HysopFFTWarning HysopFFTDataLayoutError

Methods defined here:

simd_alignment, is_byte_aligned, mk_shape, mk_view

class hysop.numerics.fft.fft.FFTI(backend, warn_on_allocation=True, error_on_allocation=False)[source]

Bases: object

Interface to compute local to process FFT-like transforms. Common inteface for all array backends, based on the numpy.fft interface.

Standard FFTs: complex to complex (C2C)

fft() Compute the 1-dimensional discrete Fourier Transform. ifft() Compute the 1-dimensional inverse discrete Fourier Transform.

Real data FFTS: real to complex (R2C) and complex to real (C2R)

rfft() Compute the 1-dimensional discrete Fourier Transform for real input. irfft() Compute the inverse of the 1-dimensional FFT of real input.

Real FFTS: real to real (R2R)

dct() Compute one of the discrete cosine transforms of a real input. dst() Compute one of the discrete sine transforms of a real input.

Supported R2R transforms are at least:

DCT-I, DCT-II, DCT-III DST-I, DST-II, DST-III

Other R2R transforms:

DCT-IV and DCT-IV are only supported by the FFTW backend at this time. DCT-V to DCT-VIII and DST-V to DST-VII are not supported by any FFT backend.

About floating point precision:

By default, both simple and double precision are supported. numpy only supports double precision (simple precision is supported by casting). FFTW also supports long double precision.

Normalization:

The default normalization has the direct transforms unscaled and the inverse transform is scaled by 1/N where N is the logical size of the transform. N should not to be confused with the physical size of the input arrays n:

FFT, RFFT: N = n DCT-II, DCT-III, DCT-IV: N = 2*n DST-II, DST-III, DST-IV: N = 2*n DCT-I: N = 2*(n-1) DST-I: N = 2*(n+1)

Orthogonal normalization is not supported by default, however the user may specify its custom normalization for each transform via the ‘scaling’ keyword parameter.

Inverse transforms (up to scaling):

Just add i in front of the method name to get the inverse transform with good scaling. For a given transform T, iT(T(X)) should always yield X within machine accuracy.

Underlying inverse transform mapping is:

FFT: IFFT RFFT: IRFFT

DCT-I: DCT-I DCT-II: DCT-III DCT-III: DCT-II DCT-IV: DCT-IV

DST-I: DST-I DST-II: DST-III DST-III: DST-II DST-IV: DST-IV

Other methods that this interface defines:

*Create queue *Transpose *Copy *Zero fill

Those methods will be used by the n-dimensional planner.

Initializes the interface and default supported real and complex types.

allocate_output(out, shape, dtype)[source]

Alocate output if required and check shape and dtype.

allocate_plans(op, plans, tmp_buffer=None)[source]

Allocate and share a buffer on given backend to a group of plans.

check_backend(backend, enable_opencl_host_buffer_mapping)[source]
abstract dct(a, out=None, type=2, axis=-1, **kwds)[source]

Compute the one-dimensional Cosine Transform of specified type.

Parameters:
  • a (array_like) – Real input array.

  • out (array_like) – Real output array of matching input type and shape.

  • axis (int, optional) – Axis over witch to compute the transform. Defaults to last axis.

Return type:

(shape, dtype) of the output array determined from the input array.

classmethod default_interface(**kwds)[source]

Get the default FFT interface.

classmethod default_interface_from_backend(backend, enable_opencl_host_buffer_mapping, **kwds)[source]
abstract dst(a, out=None, type=2, axis=-1, **kwds)[source]

Compute the one-dimensional Sine Transform of specified type.

Parameters:
  • a (array_like) – Real input array.

  • out (array_like) – Real output array of matching input type and shape.

  • axis (int, optional) – Axis over witch to compute the transform. Defaults to last axis.

Return type:

(shape, dtype) of the output array determined from the input array.

abstract fft(a, out, axis=-1, **kwds)[source]

Compute the unscaled one-dimensional complex to complex discrete Fourier Transform.

Parameters:
  • a (array_like of np.complex64 or np.complex128) – Complex input array.

  • out (array_like of np.complex64 or np.complex128) – Complex output array of the same shape and dtype as the input.

  • axis (int, optional) – Axis over witch to compute the FFT. Defaults to last axis.

Return type:

(shape, dtype) of the output array determined from the input array.

Notes

N = a.shape[axis] out[0] will contain the sum of the signal (zero-frequency term always real for real inputs).

If N is even:

out[1:N/2] contains the positive frequency terms. out[N/2] contains the Nyquist frequency (always real for real inputs). out[N/2+1:] contains the negative frequency terms.

Else if N is odd:

out[1:(N-1)/2] contains the positive frequency terms. out[(N-1)/2:] contains the negative frequency terms.

get_transform(transform)[source]
abstract idct(a, out=None, type=2, axis=-1, **kwds)[source]

Compute the one-dimensional Inverse Cosine Transform of specified type.

Default scaling is 1/(2*N) for IDCT type (2,3,4) and

1/(2*N-2) for IDCT type 1.

Parameters:
  • a (array_like) – Real input array.

  • out (array_like) – Real output array of matching input type and shape.

  • axis (int, optional) – Axis over witch to compute the transform. Defaults to last axis.

Returns:

  • (shape, dtype, inverse_type, logical_size) of the output array determined from the input

  • array.

abstract idst(a, out=None, type=2, axis=-1, **kwds)[source]

Compute the one-dimensional Inverse Sine Transform of specified type.

Default scaling is 1/(2*N) for IDST type (2,3,4) and

1/(2*N+2) for IDST type 1.

Parameters:
  • a (array_like) – Real input array.

  • out (array_like) – Real output array of matching input type and shape.

  • axis (int, optional) – Axis over witch to compute the transform. Defaults to last axis.

Returns:

  • (shape, dtype, inverse_type, logical_size) of the output array determined from the input

  • array.

abstract ifft(a, out, axis=-1, **kwds)[source]

Compute the one-dimensional complex to complex discrete Fourier Transform scaled by 1/N.

Parameters:
  • a (array_like of np.complex64 or np.complex128) – Complex input array.

  • out (array_like of np.complex64 or np.complex128) – Complex output array of the same shape and dtype as the input.

  • axis (int, optional) – Axis over witch to compute the FFT. Defaults to last axis.

Return type:

(shape, dtype, logical_size) of the output array determined from the input array.

abstract irfft(a, out, n=None, axis=-1, **kwds)[source]

Compute the one-dimensional hermitian complex to real discrete Fourier Transform scaled by 1/N.

Parameters:
  • a (array_like of np.complex64 or np.complex128) – Complex input array.

  • out (array_like of np.float32 or np.float64) –

    Real output array of matching real type. out.shape[…] = a.shape[…] Last axis should match forward transform size used:

    1. out.shape[axis] = 2*(a.shape[axis]-1)

    2. out.shape[axis] = 2*(a.shape[axis]-1) + 1

  • n (int, optional) – Length of the transformed axis of the output. ie: n should be in [2*(a.shape[axis]-1), 2*(a.shape[axis]-1)+1]

  • axis (int, optional) – Axis over witch to compute the transform. Defaults to last axis.

Notes

To get an odd number of output points, n or out must be specified.

Returns:

  • (shape, dtype, logical_size) of the output array determined from the input array,

  • out and n.

abstract new_queue(tg, name)[source]

Return a FFTQueue object valid with this backend.

abstract plan_accumulate(tg, src, dst)[source]

Plan an accumulation from src into dst.

abstract plan_compute_energy(tg, fshape, src, dst, transforms, mutexes=None)[source]

Plan to compute energy from src to energy.

abstract plan_copy(tg, src, dst)[source]

Plan a copy from src to dst.

abstract plan_fill_zeros(tg, a, slices)[source]

Plan to fill every input slices of input array a with zeroes.

abstract plan_transpose(tg, src, dst, axes)[source]

Plan a transpose from src to dst using given axes.

abstract rfft(a, out, axis=-1, **kwds)[source]

Compute the unscaled one-dimensional real to hermitian complex discrete Fourier Transform.

Parameters:
  • a (array_like of np.float32 or np.float64) – Real input array.

  • out (array_like of np.complex64 or np.complex128) – Complex output array of matching complex dtype. out.shape[…] = a.shape[…] out.shape[axis] = a.shape[axis]//2 + 1

  • axis (int, optional) – Axis over witch to compute the transform. Defaults to last axis.

Return type:

(shape, dtype) of the output array determined from the input array.

Notes

For real inputs there is no information in the negative frequency components that is not already available from the positive frequency component because of the Hermitian symmetry.

N = out.shape[axis] = a.shape[axis]//2 + 1 out[0] will contain the sum of the signal (zero-frequency term, always real). If N is even:

out[1:N/2] contains the positive frequency terms. out[N/2] contains the Nyquist frequency (always real).

Else if N is odd:

out[1:(N+1)/2] contains the positive frequency terms.

class hysop.numerics.fft.fft.FFTPlanI(verbose=False)[source]

Bases: object

Common inteface for FFT plans. Basically just a functor that holds relevant data to execute a preconfigurated FFT-like tranform.

__call__(a=None, out=None, **kwds)[source]

Apply the FFT plan to possibly new input or output arrays.

allocate(buf=None)[source]

Provide or allocate required temporary buffer.

abstract execute()[source]

Execute the FFT plan on current input and output array.

abstract input_array()[source]

Return currently planned input array.

abstract output_array()[source]

Return currently planned output array.

property required_buffer_size

Return the required temporary buffer size in bytes to compute the transform.

setup(queue=None)[source]

Method that has to be called before any call to execute.

class hysop.numerics.fft.fft.FFTQueueI[source]

Bases: object

Command queue like objects to define n-dimensional transforms.

__call__(**kwds)[source]

Alias for execute.

abstract __iadd__(*plans)[source]

Add a plan to the queue.

abstract execute(wait_for=None)[source]

Execute all planned plans.

exception hysop.numerics.fft.fft.HysopFFTDataLayoutError[source]

Bases: ValueError

Specific error to raise for incompatible strides.

exception hysop.numerics.fft.fft.HysopFFTWarning[source]

Bases: HysopWarning

Specific tag to issue FFT related warnings.

hysop.numerics.fft.fft.is_byte_aligned(array)[source]
hysop.numerics.fft.fft.mk_shape(base_shape, axis, N)[source]

Utility function to create a modified shape from base_shape on a specified axis.

hysop.numerics.fft.fft.mk_view(ndim, axis, *args, **kwds)[source]

Utility function to create a view on a n-dimensional array. Returns a tuple containing default axe view (ellipsis) on all axis but the one specified by ‘axis’ which is replaced by slice(*args) if len(args)>1 else args[0].